GeoLib Module

History

current version 0.9 - 16 July 2018

Version Date Comment
0.1 22/Dec/2008 Original code. giovanni ravazzani
0.2 09/Aug/2009 Added ScanDatum to convert string to numeric code
0.3 12/Jul/2010 Support for Hotine Oblique Mercator
0.4 15/Jul/2010 Support for Swiss Oblique Cylindrical
0.5 25/Nov/2010 Pass to double precision in ConvertTransverseMercatorToGeodetic
0.6 03/Feb/2011 Added function for computing distance between two points
0.7 20/May/2012 Extend Coordinate type with elevation element
0.8 13/Jan/2016 DecodeEPSG to read epsg code and set CRS
0.9 16/Jul/2018 Added functions to compute the bearing and direction angle between two points

License

license: GNU GPL http://www.gnu.org/licenses/

This file is part of

MOSAICO -- MOdular library for raSter bAsed hydrologIcal appliCatiOn.

Copyright (C) 2011 Giovanni Ravazzani

Module Description

set of fortran routines that allows you to convert geographic coordinates among a variety of coordinate systems, map projections, and datums.
The function Convert is the heart of the GeoLib module. It is this function that actually performs all coordinate conversions and datum shifts. Given the input and output datums, the input and output coordinate reference system, including the values of any projection parameters, and the input coordinates, it produces the output coordinates. Methods: The general case is handled in three stages: 1. Convert the input coordinates from the input coordinate reference system to geodetic, 2. Shift the intermediate geodetic coordinates from the input datum to the output datum, converting between ellipsoid if necessary, 3. Convert the shifted intermediate geodetic coordinates to the output coordinate reference system. The first stage is accomplished by a case statement on the input coordinate reference system, with a case for each of the coordinate reference systems supported. If the input coordinates are already in geodetic, they are simply copied to the intermediate coordinates. The second stage is accomplished in two steps. First, the intermediate Geodetic coordinates are shifted to WGS 84, if necessary. Second, the shifted intermediate geodetic coordinates are shifted to the output datum. The third stage is similar to the first stage. A case statement on the output coordinate reference system, with a case for each of the coordinate reference frames supported, is used to convert the shifted intermediate Geodetic coordinates to the output coordinate reference system. If the output coordinate reference system is Geodetic, the shifted intermediate coordinates are simply copied to the output coordinates.


References:

GeoTrans: http://earth-info.nga.mil/GandG/geotrans/index.html

European Petroleum Survey Group (EPSG): http://www.epsg.org

http://spatialreference.org/

http://www.borneo.name/topo/roma40.html



Variables

Type Visibility Attributes Name Initial
integer, public, parameter :: CH1903 = 6149

World Geodetic System 1984 European Datum 1950 Monte Mario Italy Switzerland

integer, public, parameter :: EAST = 3
integer, public, parameter :: ED50 = 6230

World Geodetic System 1984 European Datum 1950 Monte Mario Italy Switzerland

integer, public, parameter :: GAUSS_BOAGA = 3
integer, public, parameter :: GEODETIC = 1
integer, public, parameter :: HOM = 5
integer, public, parameter :: NORTH = 1
integer, public, parameter :: ROME40 = 6265

World Geodetic System 1984 European Datum 1950 Monte Mario Italy Switzerland

integer, public, parameter :: SOC = 6
integer, public, parameter :: SOUTH = 2
integer, public, parameter :: TM = 4
integer, public, parameter :: UTM = 2
integer, public, parameter :: WEST = 4
integer, public, parameter :: WGS84 = 6326

World Geodetic System 1984 European Datum 1950 Monte Mario Italy Switzerland

type(Coordinate), public :: point1

can be used to calculate distance between two points

type(Coordinate), public :: point2

can be used to calculate distance between two points

integer, private, parameter :: GB_centM = 2

Zone override flag Zone override flag

integer, private, parameter :: GB_false_easting = 4

Zone override flag Zone override flag

integer, private, parameter :: GB_false_northing = 5

Zone override flag Zone override flag

integer, private, parameter :: GB_fuse = 3

Zone override flag Zone override flag

integer, private, parameter :: GB_lat0 = 1

Zone override flag Zone override flag

integer, private, parameter :: GB_override = 7

Zone override flag Zone override flag

integer, private, parameter :: GB_scale_factor = 6

Zone override flag Zone override flag

integer, private, parameter :: GEO_a = 1

Zone override flag Zone override flag

integer, private, parameter :: GEO_b = 2

Zone override flag Zone override flag

integer, private, parameter :: GEO_invf = 3

Zone override flag Zone override flag

integer, private, parameter :: GEO_prime = 4

Zone override flag Zone override flag

integer, private, parameter :: HOM_azimuth = 3

Zone override flag Zone override flag

integer, private, parameter :: HOM_centM = 2

Zone override flag Zone override flag

integer, private, parameter :: HOM_false_easting = 4

Zone override flag Zone override flag

integer, private, parameter :: HOM_false_northing = 5

Zone override flag Zone override flag

integer, private, parameter :: HOM_lat0 = 1

Zone override flag Zone override flag

integer, private, parameter :: HOM_scale_factor = 6

Zone override flag Zone override flag

integer, private, parameter :: SOC_azimuth = 3

Zone override flag Zone override flag

integer, private, parameter :: SOC_false_easting = 4

Zone override flag Zone override flag

integer, private, parameter :: SOC_false_northing = 5

Zone override flag Zone override flag

integer, private, parameter :: SOC_latc = 1

Zone override flag Zone override flag

integer, private, parameter :: SOC_lonc = 2

Zone override flag Zone override flag

integer, private, parameter :: SOC_scale_factor = 6

Zone override flag Zone override flag

integer, private, parameter :: TM_centM = 2

Zone override flag Zone override flag

integer, private, parameter :: TM_false_easting = 3

Zone override flag Zone override flag

integer, private, parameter :: TM_false_northing = 4

Zone override flag Zone override flag

integer, private, parameter :: TM_lat0 = 1

Zone override flag Zone override flag

integer, private, parameter :: TM_scale_factor = 5

Zone override flag Zone override flag

integer, private, parameter :: UTM_centM = 2

Zone override flag Zone override flag

integer, private, parameter :: UTM_false_easting = 5

Zone override flag Zone override flag

integer, private, parameter :: UTM_false_northing = 6

Zone override flag Zone override flag

integer, private, parameter :: UTM_hemisphere = 4

Zone override flag Zone override flag

integer, private, parameter :: UTM_lat0 = 1

Zone override flag Zone override flag

integer, private, parameter :: UTM_override = 8

Zone override flag Zone override flag

integer, private, parameter :: UTM_scale_factor = 7

Zone override flag Zone override flag

integer, private, parameter :: UTM_zone = 3

Zone override flag Zone override flag

type(datum), private, ALLOCATABLE :: datums(:)

array of available datums

type(Ellipsoid), private, ALLOCATABLE :: ellps(:)

array of available ellipsoids

real(kind=Float), private, parameter :: null_float = -9999.9
character(len=0), private, parameter :: null_string = ''

Interfaces

public interface ASSIGNMENT( = )

  • private subroutine CopyEllipsoid(ell2, ell1)

    Create an exact copy of an ellipsoid.

    Arguments

    Type IntentOptional Attributes Name
    type(Ellipsoid), intent(out) :: ell2
    type(Ellipsoid), intent(in) :: ell1
  • private subroutine SearchEllipsoidByCode(ell, code)

    Search for ellipsoid parameters using a code.

    Arguments

    Type IntentOptional Attributes Name
    type(Ellipsoid), intent(out) :: ell
    integer, intent(in) :: code
  • private subroutine CopyDatum(datum2, datum1)

    Create an exact copy of a datum.

    Arguments

    Type IntentOptional Attributes Name
    type(datum), intent(out) :: datum2
    type(datum), intent(in) :: datum1
  • private subroutine SearchDatumByCode(dat, code)

    Search for datum parameters using a code.

    Arguments

    Type IntentOptional Attributes Name
    type(datum), intent(out) :: dat
    integer, intent(in) :: code
  • private subroutine CopyCRS(CRSout, CRSin)

    copy the content of a CRS variable in another CRS variable

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(out) :: CRSout
    type(CRS), intent(in) :: CRSin

public interface OPERATOR (==)

  • private function CRSisEqual(CRS1, CRS2) result(isEqual)

    return .TRUE. if the two coordinate reference systems are equal

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(in) :: CRS1
    type(CRS), intent(in) :: CRS2

    Return Value logical

  • private function EllipsoidIsEqual(ellps1, ellps2) result(isEqual)

    return .TRUE. if the two ellipsoids are equal

    Arguments

    Type IntentOptional Attributes Name
    type(Ellipsoid), intent(in) :: ellps1
    type(Ellipsoid), intent(in) :: ellps2

    Return Value logical

  • private function DatumIsEqual(datum1, datum2) result(isEqual)

    return .TRUE. if the two datums are equal

    Arguments

    Type IntentOptional Attributes Name
    type(datum), intent(in) :: datum1
    type(datum), intent(in) :: datum2

    Return Value logical

public interface SetCRS

  • private subroutine SetCRScoord(CRStype, datumType, coord)

    Initialize Coordinate Reference System, allocate memory and set parameters to default value if necessary. subroutine receives as input a Coordinate type argument

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: CRStype
    integer, intent(in) :: datumType
    type(Coordinate), intent(inout) :: coord
  • private subroutine SetCRSsystem(CRStype, datumType, rs)

    Initialize Coordinate Reference System, allocate memory and set parameters to default value if necessary. Subroutine receives as input a CRS type argument

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: CRStype
    integer, intent(in) :: datumType
    type(CRS), intent(inout) :: rs

    reference system

public interface SetGaussBoagaParameters

  • private subroutine SetGaussBoagaparametersSystem(system, fuse, override)

    Set parameters for Gauss Boaga reference system

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    integer(kind=short), intent(in) :: fuse
    integer(kind=short), intent(in), optional :: override
  • private subroutine SetGaussBoagaparametersCoord(coord, fuse, override)

    Set parameters for Gauss Boaga reference system

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    integer(kind=short), intent(in) :: fuse
    integer(kind=short), intent(in), optional :: override

public interface SetGeodeticParameters

  • private subroutine SetGeodeticParametersSystem(system, prime_meridian)

    Set parameters for Geodetic reference system

    longitude, with respect to Greenwich, of the prime meridian

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    real(kind=float), intent(in) :: prime_meridian
  • private subroutine SetGeodeticParametersCoord(coord, prime_meridian)

    Set parameters for Geodetic reference system

    longitude, with respect to Greenwich, of the prime meridian

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    real(kind=float), intent(in) :: prime_meridian

public interface SetHotineParameters

  • private subroutine SetHotineParametersSystem(system, lat0, centM, azimuth, falseE, falseN, k)

    Set parameters for Hotine Oblique Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    real(kind=float), intent(in) :: lat0

    Latitude of projection center [rad]

    real(kind=float), intent(in) :: centM

    Longitude of projection center [rad]

    real(kind=float), intent(in) :: azimuth

    azimuth of centerline

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

  • private subroutine SetHotineParametersCoord(coord, lat0, centM, azimuth, falseE, falseN, k)

    Set parameters for Hotine Oblique Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    real(kind=float), intent(in) :: lat0

    Latitude of projection center [rad]

    real(kind=float), intent(in) :: centM

    Longitude of projection center [rad]

    real(kind=float), intent(in) :: azimuth

    azimuth of centerline

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

public interface SetSwissParameters

  • private subroutine SetSwissParametersSystem(system, latc, lonc, azimuth, falseE, falseN, k)

    Set parameters for Swiss Oblique Cylindrical reference system

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    real(kind=float), intent(in) :: latc

    Latitude of projection center [rad]

    real(kind=float), intent(in) :: lonc

    Longitude of projection center [rad]

    real(kind=float), intent(in) :: azimuth

    azimuth of centerline

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

  • private subroutine SetSwissParametersCoord(coord, latc, lonc, azimuth, falseE, falseN, k)

    Set parameters for Swiss Oblique Cylindrical reference system

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    real(kind=float), intent(in) :: latc

    Latitude of projection center [rad]

    real(kind=float), intent(in) :: lonc

    Longitude of projection center [rad]

    real(kind=float), intent(in) :: azimuth

    azimuth of centerline

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

  • private subroutine SetTransverseMercatorParametersSystem(system, lat0, centM, falseE, falseN, k)

    Set parameters for Transverse Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    real(kind=float), intent(in) :: lat0

    Latitude in radians at the origin of the projection

    real(kind=float), intent(in) :: centM

    Longitude in radians at the center of the projection

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

  • private subroutine SetTransverseMercatorParametersCoord(coord, lat0, centM, falseE, falseN, k)

    Set parameters for Transverse Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    real(kind=float), intent(in) :: lat0

    Latitude in radians at the origin of the projection

    real(kind=float), intent(in) :: centM

    Longitude in radians at the center of the projection

    real(kind=float), intent(in) :: falseE

    Easting/X at the center of the projection

    real(kind=float), intent(in) :: falseN

    Northing/Y at the center of the projection

    real(kind=float), intent(in) :: k

    scale factor

public interface SetUTMparameters

  • private subroutine SetUTMparametersSystem(system, zone, hemisphere, override)

    Set parameters for Universal Transverse Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(CRS), intent(inout) :: system
    integer(kind=short), intent(in) :: zone
    integer(kind=short), intent(in) :: hemisphere
    integer(kind=short), intent(in), optional :: override
  • private subroutine SetUTMparametersCoord(coord, zone, hemisphere, override)

    Set parameters for Universal Transverse Mercator reference system

    Arguments

    Type IntentOptional Attributes Name
    type(Coordinate), intent(inout) :: coord
    integer(kind=short), intent(in) :: zone
    integer(kind=short), intent(in) :: hemisphere
    integer(kind=short), intent(in), optional :: override

Derived Types

type, public ::  CRS

Components

Type Visibility Attributes Name Initial
integer(kind=short), public :: basic

number of basic parameters for the reference system

type(datum), public :: datum
character(len=100), public, ALLOCATABLE :: description(:)

description of the parameter (eg. zone for UTM, etc..)

type(Ellipsoid), public :: ellipsoid
integer(kind=short), public :: epsg

epsg code

character(len=50), public :: name

name of CRS

real(kind=FLOAT), public, ALLOCATABLE :: param(:)

required parameters for the definition of reference system

integer(kind=short), public :: system

geodetic, UTM, Gaus Boaga, etc..

type, public ::  Coordinate

Components

Type Visibility Attributes Name Initial
real(kind=Float), public :: easting

X coordinate, longitude

real(kind=Float), public :: elevation

Z coordinate

real(kind=Float), public :: northing

Y coordinate, latitude

type(CRS), public :: system

coordinate reference system

type, public ::  Ellipsoid

Components

Type Visibility Attributes Name Initial
real(kind=Float), public :: a

semi-major axis

real(kind=Float), public :: b

semi-minor axis

integer(kind=short), public :: code

EPSG code

real(kind=Float), public :: e

eccentricity

real(kind=Float), public :: e_second

second eccentricity

character(len=100), public :: epsg

EPSG string

real(kind=Float), public :: f

flattening

real(kind=Float), public :: inv_f

a/(a-b)

character(len=100), public :: name

type, public ::  datum

Components

Type Visibility Attributes Name Initial
integer(kind=short), public :: code
real(kind=Float), public :: dx

entity of x-axis translation

real(kind=Float), public :: dy

entity of y-axis translation

real(kind=Float), public :: dz

entity of z-axis translation

integer(kind=short), public :: ellipsoid

reference ellipsoid code

character(len=100), public :: epsg

EPSG string

integer(kind=short), public :: method

EPSG code defining operation method for translation to WGS84

character(len=100), public :: name
integer(kind=short), public :: prime_meridian

EPSG code defining prime meridian


Functions

public function BearingAngle(pointA, pointB) result(angle)

Find the bearing angle (radians) between two points in a 2D space, as the angle measured in the clockwise direction from the north line with A as the origin to the line segment AB. Assume radians as input unit for geodetic reference system.

Read more…

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: pointA

origin point

type(Coordinate), intent(in) :: pointB

ending point

Return Value real(kind=float)

bearing angle (radians)

public function DecodeEPSG(epsg) result(gridmapping)

get EPSG code as input and set coordinate reference system

Arguments

Type IntentOptional Attributes Name
integer(kind=short), intent(in) :: epsg

Return Value type(CRS)

public function DirectionAngle(pointA, pointB) result(angle)

Find the direction angle (radians) between two points in a 2D space, measured in the clockwise direction. Formula gives values 0°<θ<90° for lines with positive slope and values 90°<θ<180° for lines with negative slope. The result does depend on which point is point 1 and which is point 2. Assume radians as input unit for geodetic reference system.

Read more…

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: pointA

point A

type(Coordinate), intent(in) :: pointB

point B

Return Value real(kind=float)

bearing angle (radians)

public function Distance(point1, point2) result(distance)

compute distance between two points in, either, projected or geodetic coordinate reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: point1
type(Coordinate), intent(in) :: point2

Return Value real(kind=double)

public function ScanDatum(datumString) result(datumCode)

return datum numeric code given a string

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: datumString

Return Value integer

private function CRSisEqual(CRS1, CRS2) result(isEqual)

return .TRUE. if the two coordinate reference systems are equal

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(in) :: CRS1
type(CRS), intent(in) :: CRS2

Return Value logical

private function DatumIsEqual(datum1, datum2) result(isEqual)

return .TRUE. if the two datums are equal

Arguments

Type IntentOptional Attributes Name
type(datum), intent(in) :: datum1
type(datum), intent(in) :: datum2

Return Value logical

private function EllipsoidIsEqual(ellps1, ellps2) result(isEqual)

return .TRUE. if the two ellipsoids are equal

Arguments

Type IntentOptional Attributes Name
type(Ellipsoid), intent(in) :: ellps1
type(Ellipsoid), intent(in) :: ellps2

Return Value logical

private function MeridionalDistance(lat, a, e) result(M)

Compute the distance in meters along a meridian from equator to given latitude.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: lat

latitude in [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

exentricity

Return Value real(kind=float)

meridional distance [m]

private function Nu(lat, a, e) result(radius)

Compute the radius of curvature of the ellipsoid perpendicular to the meridian at given latitude

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: lat

latitude [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

Return Value real(kind=float)

private function Rho(lat, a, e) result(radius)

Compute the radius of curvature of the ellipsoid in the plane of the meridian at given latitude

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: lat

latitude [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

exentricity

Return Value real(kind=float)


Subroutines

public subroutine Convert(input, output)

The subroutine converts the current input coordinates in the coordinate system defined by the current input coordinate system parameters and input datum, into output coordinates in the coordinate system defined by the output coordinate system parameters and output datum.

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: input
type(Coordinate), intent(inout) :: output

public subroutine ConvertGeodeticToSwiss(lon, lat, k, lon0, lat0, azimuth, a, e, falseN, falseE, x, y)

The subroutine converts geodetic (latitude and longitude) to Swiss Oblique Cylindrical projection(easting and northing) coordinates

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: lon

geodetic longitude [radians]

real(kind=float), intent(in) :: lat

geodetic latitude [radians]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: lon0

longitude of center [radians]

real(kind=float), intent(in) :: lat0

latitude of center [radians]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: x

easting coordinate [m]

real(kind=float), intent(out) :: y

northing coordinate [m]

public subroutine GeoInit(file)

Initialize parameters for conversion

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: file

public subroutine SetCoord(x, y, coord)

assign easting and northing coordinates

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate

real(kind=float), intent(in) :: y

northing coordinate

type(Coordinate), intent(inout) :: coord

private subroutine ConvertGaussBoagaToGeodetic(x, y, k, centM, lat0, a, e, eb, falseN, falseE, lon, lat)

The subroutine converts Gauss Boaga projection for Italy (easting and northing) coordinates to geodetic (latitude and longitude) coordinates

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate [m]

real(kind=float), intent(in) :: y

northing coordinate [m]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: lon

geodetic longitude [radians]

real(kind=float), intent(out) :: lat

geodetic latitude [radians]

private subroutine ConvertGeodeticToGaussBoaga(lon, lat, k, centM, lat0, a, e, eb, falseN, falseE, x, y)

The subroutine converts geodetic(latitude and longitude) coordinates to Gauss Boaga (easting and northing) coordinates, according to ROME40 datum for Italy (MonteMario) coordinates.

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(inout) :: lon

geodetic longitude [radians]

real(kind=float), intent(inout) :: lat

geodetic latitude [radians]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: x

easting coordinate [m]

real(kind=float), intent(out) :: y

northing coordinate [m]

private subroutine ConvertGeodeticToHotine(lon, lat, k, lon0, lat0, azimuth, a, e, falseN, falseE, x, y)

The subroutine converts geodetic (latitude and longitude) coordinates to Hotine Oblique Mercator projection (easting and northing)

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: lon

geodetic longitude [radians]

real(kind=float), intent(in) :: lat

geodetic latitude [radians]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: lon0

longitude of center [radians]

real(kind=float), intent(in) :: lat0

latitude of center [radians]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: x

easting coordinate [m]

real(kind=float), intent(out) :: y

northing coordinate [m]

private subroutine ConvertGeodeticToTransverseMercator(lon, lat, k, centM, lat0, a, e, eb, falseN, falseE, x, y)

The subroutine converts geodetic(latitude and longitude) coordinates to Transverse Mercator projection (easting and northing) coordinates, according to the current ellipsoid and Transverse Mercator projection coordinates.

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(inout) :: lon

geodetic longitude [radians]

real(kind=float), intent(in) :: lat

geodetic latitude [radians]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: x

easting coordinate [m]

real(kind=float), intent(out) :: y

northing coordinate [m]

private subroutine ConvertGeodeticToUTM(lon, lat, k, centM, lat0, a, e, eb, falseN, falseE, x, y, zone, hemisphere, override)

The subroutine converts geodetic(latitude and longitude) coordinates to UTM projection (easting and northing) coordinates, according to the current ellipsoid and UTM projection coordinates.

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(inout) :: lon

geodetic longitude [radians]

real(kind=float), intent(inout) :: lat

geodetic latitude [radians]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(inout) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(inout) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: x

easting coordinate [m]

real(kind=float), intent(out) :: y

northing coordinate [m]

real(kind=float), intent(inout) :: zone

zone recalculated if override = 0.

real(kind=float), intent(inout) :: hemisphere

hemisphere recalculated if override = 0.

real(kind=float), intent(in) :: override

override option

private subroutine ConvertHotineToGeodetic(x, y, k, lon0, lat0, azimuth, a, e, falseN, falseE, lon, lat)

The subroutine converts Hotine Oblique Mercator projection (easting and northing) coordinates to geodetic (latitude and longitude) coordinates

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate [m]

real(kind=float), intent(in) :: y

northing coordinate [m]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: lon0

longitude of center [radians]

real(kind=float), intent(in) :: lat0

latitude of center [radians]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: lon

geodetic longitude [radians]

real(kind=float), intent(out) :: lat

geodetic latitude [radians]

private subroutine ConvertSwissToGeodetic(x, y, k, lon0, lat0, azimuth, a, e, falseN, falseE, lon, lat)

The subroutine converts Swiss Oblique Cylindrical projection (easting and northing) coordinates to geodetic (latitude and longitude) coordinates

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate [m]

real(kind=float), intent(in) :: y

northing coordinate [m]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: lon0

longitude of center [radians]

real(kind=float), intent(in) :: lat0

latitude of center [radians]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: lon

geodetic longitude [radians]

real(kind=float), intent(out) :: lat

geodetic latitude [radians]

private subroutine ConvertTransverseMercatorToGeodetic(x, y, k, centM, lat0, a, e, eb, falseN, falseE, lon, lat)

The subroutine converts Transverse Mercator projection (easting and northing) coordinates to geodetic (latitude and longitude) coordinates, according to the current ellipsoid and Transverse Mercator projection parameters.

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate [m]

real(kind=float), intent(in) :: y

northing coordinate [m]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(out) :: lon

geodetic longitude [radians]

real(kind=float), intent(out) :: lat

geodetic latitude [radians]

private subroutine ConvertUTMtoGeodetic(x, y, k, centM, lat0, a, e, eb, falseN, falseE, override, lon, lat)

The subroutine converts Universal Transverse Mercator projection (easting and northing) coordinates to geodetic (latitude and longitude) coordinates, according to the current ellipsoid and UTM projection parameters.

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: x

easting coordinate [m]

real(kind=float), intent(in) :: y

northing coordinate [m]

real(kind=float), intent(in) :: k

scale factor

real(kind=float), intent(in) :: centM

central meridian [radians]

real(kind=float), intent(in) :: lat0

latitude of origin [radians]

real(kind=float), intent(in) :: a

semimajor axis [m]

real(kind=float), intent(in) :: e

eccentricity

real(kind=float), intent(in) :: eb

second eccentricity

real(kind=float), intent(in) :: falseN

false northing

real(kind=float), intent(in) :: falseE

false easting

real(kind=float), intent(in) :: override

override option

real(kind=float), intent(out) :: lon

geodetic longitude [radians]

real(kind=float), intent(out) :: lat

geodetic latitude [radians]

private subroutine CopyCRS(CRSout, CRSin)

copy the content of a CRS variable in another CRS variable

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(out) :: CRSout
type(CRS), intent(in) :: CRSin

private subroutine CopyDatum(datum2, datum1)

Create an exact copy of a datum.

Arguments

Type IntentOptional Attributes Name
type(datum), intent(out) :: datum2
type(datum), intent(in) :: datum1

private subroutine CopyEllipsoid(ell2, ell1)

Create an exact copy of an ellipsoid.

Arguments

Type IntentOptional Attributes Name
type(Ellipsoid), intent(out) :: ell2
type(Ellipsoid), intent(in) :: ell1

private subroutine GeodeticShiftFromWGS84(input, WGS84lat, WGS84lon, latOut, lonOut)

shifts geodetic coordinates relative to WGS84 to geodetic coordinates relative to a given local datum.

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: input
real(kind=float), intent(in) :: WGS84lat
real(kind=float), intent(in) :: WGS84lon
real(kind=float), intent(out) :: latOut
real(kind=float), intent(out) :: lonOut

private subroutine GeodeticShiftToWGS84(input, latIn, lonIn, WGS84lat, WGS84lon)

shifts geodetic coordinates relative to a given source datum to geodetic coordinates relative to WGS84.

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(in) :: input
real(kind=float), intent(in) :: latIn
real(kind=float), intent(in) :: lonIn
real(kind=float), intent(out) :: WGS84lat
real(kind=float), intent(out) :: WGS84lon

private subroutine MolodenskyShift(a, da, f, df, dx, dy, dz, Lat_in, Lon_in, Hgt_in, WGS84_Lat, WGS84_Lon, WGS84_Hgt)

shifts geodetic coordinates using the Molodensky method

Arguments

Type IntentOptional Attributes Name
real(kind=float), intent(in) :: a
real(kind=float), intent(in) :: da
real(kind=float), intent(in) :: f
real(kind=float), intent(in) :: df
real(kind=float), intent(in) :: dx
real(kind=float), intent(in) :: dy
real(kind=float), intent(in) :: dz
real(kind=float), intent(in) :: Lat_in
real(kind=float), intent(in) :: Lon_in
real(kind=float), intent(in) :: Hgt_in
real(kind=float), intent(out) :: WGS84_Lat
real(kind=float), intent(out) :: WGS84_Lon
real(kind=float), intent(out) :: WGS84_Hgt

private subroutine SearchDatumByCode(dat, code)

Search for datum parameters using a code.

Arguments

Type IntentOptional Attributes Name
type(datum), intent(out) :: dat
integer, intent(in) :: code

private subroutine SearchEllipsoidByCode(ell, code)

Search for ellipsoid parameters using a code.

Arguments

Type IntentOptional Attributes Name
type(Ellipsoid), intent(out) :: ell
integer, intent(in) :: code

private subroutine SetCRScoord(CRStype, datumType, coord)

Initialize Coordinate Reference System, allocate memory and set parameters to default value if necessary. subroutine receives as input a Coordinate type argument

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: CRStype
integer, intent(in) :: datumType
type(Coordinate), intent(inout) :: coord

private subroutine SetCRSsystem(CRStype, datumType, rs)

Initialize Coordinate Reference System, allocate memory and set parameters to default value if necessary. Subroutine receives as input a CRS type argument

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: CRStype
integer, intent(in) :: datumType
type(CRS), intent(inout) :: rs

reference system

private subroutine SetGaussBoagaparametersCoord(coord, fuse, override)

Set parameters for Gauss Boaga reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
integer(kind=short), intent(in) :: fuse
integer(kind=short), intent(in), optional :: override

private subroutine SetGaussBoagaparametersSystem(system, fuse, override)

Set parameters for Gauss Boaga reference system

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
integer(kind=short), intent(in) :: fuse
integer(kind=short), intent(in), optional :: override

private subroutine SetGeodeticParametersCoord(coord, prime_meridian)

Set parameters for Geodetic reference system

Read more…

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
real(kind=float), intent(in) :: prime_meridian

private subroutine SetGeodeticParametersSystem(system, prime_meridian)

Set parameters for Geodetic reference system

Read more…

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
real(kind=float), intent(in) :: prime_meridian

private subroutine SetHotineParametersCoord(coord, lat0, centM, azimuth, falseE, falseN, k)

Set parameters for Hotine Oblique Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
real(kind=float), intent(in) :: lat0

Latitude of projection center [rad]

real(kind=float), intent(in) :: centM

Longitude of projection center [rad]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetHotineParametersSystem(system, lat0, centM, azimuth, falseE, falseN, k)

Set parameters for Hotine Oblique Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
real(kind=float), intent(in) :: lat0

Latitude of projection center [rad]

real(kind=float), intent(in) :: centM

Longitude of projection center [rad]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetSwissParametersCoord(coord, latc, lonc, azimuth, falseE, falseN, k)

Set parameters for Swiss Oblique Cylindrical reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
real(kind=float), intent(in) :: latc

Latitude of projection center [rad]

real(kind=float), intent(in) :: lonc

Longitude of projection center [rad]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetSwissParametersSystem(system, latc, lonc, azimuth, falseE, falseN, k)

Set parameters for Swiss Oblique Cylindrical reference system

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
real(kind=float), intent(in) :: latc

Latitude of projection center [rad]

real(kind=float), intent(in) :: lonc

Longitude of projection center [rad]

real(kind=float), intent(in) :: azimuth

azimuth of centerline

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetTransverseMercatorParametersCoord(coord, lat0, centM, falseE, falseN, k)

Set parameters for Transverse Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
real(kind=float), intent(in) :: lat0

Latitude in radians at the origin of the projection

real(kind=float), intent(in) :: centM

Longitude in radians at the center of the projection

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetTransverseMercatorParametersSystem(system, lat0, centM, falseE, falseN, k)

Set parameters for Transverse Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
real(kind=float), intent(in) :: lat0

Latitude in radians at the origin of the projection

real(kind=float), intent(in) :: centM

Longitude in radians at the center of the projection

real(kind=float), intent(in) :: falseE

Easting/X at the center of the projection

real(kind=float), intent(in) :: falseN

Northing/Y at the center of the projection

real(kind=float), intent(in) :: k

scale factor

private subroutine SetUTMparametersCoord(coord, zone, hemisphere, override)

Set parameters for Universal Transverse Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(Coordinate), intent(inout) :: coord
integer(kind=short), intent(in) :: zone
integer(kind=short), intent(in) :: hemisphere
integer(kind=short), intent(in), optional :: override

private subroutine SetUTMparametersSystem(system, zone, hemisphere, override)

Set parameters for Universal Transverse Mercator reference system

Arguments

Type IntentOptional Attributes Name
type(CRS), intent(inout) :: system
integer(kind=short), intent(in) :: zone
integer(kind=short), intent(in) :: hemisphere
integer(kind=short), intent(in), optional :: override